BigDFT.InputActions module
Actions to define on the Input parameters.
This module defines some of the most common actions that a BigDFT user might
like to perform on the input file. Such module therefore sets some of the keys
of the input dictionary to the values needed to perform the operations.
Users might also inspire to the actions performed in order to customize the
runs in a different way. All the functions of this module have as first
argument inp
, the dictionary of the input parameters.
Many other actions are available in BigDFT code. This module only regroups the
most common. Any of these functionalities might be removed from the input file
by the remove()
function.
Note
Any of the action of this module, including the remove()
function,
can be also applied to an instance of the
BigDFT.Inputfiles.Inputfile
class, by removing the first
argument (inp
). This adds extra flexibility as the same method may be
used to a dictionary instance or to a BigDFT input files.
See the example Example.
We now list the available methods, in order of category.
Basic Setup and common functionalities
|
Set the exchange and correlation approximation |
|
Set the wavelet grid spacing. |
|
Set the wavelet grid extension by modifying the multiplicative radii. |
|
Constrain the number of grid points in each direction. |
|
Optimize the geometry of the system |
|
Add a collinear spin polarization to the system. |
|
Charge the system |
|
Charge the system by removing one electron. Assume that the original |
|
Set the symmetry detection for the charge density and the ionic forces and stressdef set_symmetry(inp,yes=True): |
|
Apply an external electric field on the system |
|
Insert |
|
Extract a given number of empty states after the scf cycle. |
|
Define the electronic temperature, in AU ( |
|
Add an implicit solvent around the system with the |
Add Grimme's D3 correction to the energy and the forces |
Self-Consistent-Field setup and algorithms
|
Set the tolerance acceptance level for stopping the self-consistent iterations. |
Input orbitals are initialized as random coefficients |
|
|
Activates the linear scaling mode |
|
Set the algorithm for scf. |
|
Method for guessing the kernel at restart. |
|
In the linear scaling mode, use NTPoly as the solver. |
|
Methods for scf of the density kernel. |
|
Methods for scf of the kernel coefficients. |
|
Methods for scf of the Support functions. |
Input-Output and restart
|
Set the code to write the orbitals on disk in the provided format |
|
Read the orbitals from data directory, if available. |
Write the charge density on the disk after the last scf convergence |
|
|
Modify the name of the |
|
Associate the data of the run of a given logfile to the input by retrieving the data directory name of the logfile. |
|
Write the matrices of the linear scaling formats. |
|
Writes the |
|
Dump the support functions. |
Post-Processing and other functionalities
|
Extract the dipole moment from the total charge density. |
|
Calculate the Partial Density of states. |
|
Perform a Casida TDDFT coupling matrix extraction. |
|
Calculate the multipoles on the atoms based on SF. |
|
Include a cdft constraint into the SCF cycle. |
Setting atomic-based information and other functionalities
|
Insert the atomic positions as a part of the input dictionary |
|
Employ gpu acceleration when available, for convolutions and Fock operator |
|
Control the occupation number of the KS orbitals. |
|
Set the external potential to which the system is submitted outside the QM region |
|
Set the K-point mesh. |
|
Employs all the "psppar.*" files of the given directory as pseudopotentials |
|
Employ the given PSP file for the provided element |
|
Employed the built in Non Linear Core Correction (NLCC) pseudopotentials. |
|
Load a profile or list of profiles. |
|
Remove action from the input dictionary. |
Note
Each of the actions here must have default value for the arguments
(except the input dictionary inp
). This is needed for a good behaviour
of the function remove.
We list here the extended documentation in alphabetic order.
- remove(inp, action)[source]
Remove action from the input dictionary.
Remove an action from the input file, thereby restoring the default value, as if the action were not specified.
- Parameters:
inp (dict) – dictionary to remove the action from.
action (func) – one of the actions of this module. It does not need to be
before (specified) –
effect. (in which case it produces no) –
Example
>>> from Calculators import SystemCalculator as C >>> code=C() >>> inp={} >>> set_xc(inp,'PBE') >>> write_orbitals_on_disk(inp) >>> log=code.run(input=inp) # perform calculations >>> remove(inp, write_orbitals_on_disk) #remove the action >>> read_orbitals_from_disk(inp) >>> # this will restart the scf from the previous orbitals >>> log2=code.run(input=inp)
- set_scf_convergence(inp, gnrm='default', rpnrm='default')[source]
Set the tolerance acceptance level for stopping the self-consistent iterations. Useful both for LS and CS
- set_rmult(inp, rmult=None, coarse=5.0, fine=8.0)[source]
Set the wavelet grid extension by modifying the multiplicative radii.
- Parameters:
rmult (float,list) – list of two values that have to be used for the coarse and the fine resolution grid. It may also be a scalar.
coarse (float) – if the argument
rmult
is not provided it sets the coarse radius multiplierfine (float) – if the argument
rmult
is not provided it sets the fine radius multiplier
- set_symmetry(inp, yes=True)[source]
Set the symmetry detection for the charge density and the ionic forces and stressdef set_symmetry(inp,yes=True):
- Parameters:
yes (bool) – If
False
the symmetry detection is disabled
- set_ntpoly(inp, thresh_dens=1e-06, conv_dens=0.0001, thresh_ovlp=1e-07, conv_ovlp=0.0001)[source]
In the linear scaling mode, use NTPoly as the solver.
- Parameters:
thresh_dens (float) – the threshold for filtering sparse matrices
density. (when solving for the) –
thresh_ovlp (float) – the threshold for filtering sparse matrices
inverse. (when solving for the overlap) –
conv_dens (float) – the energy value to consider the density converged.
conv_ovlp (float) – the value (in terms of norm) to consider the
converged. (overlap matrix) –
solver (int) – which solver to use; (1) TRS4, (2) TRS2,
eigensolver ((3) Dense) –
- set_mesh_sizes(inp, ngrids=64)[source]
Constrain the number of grid points in each direction. This is useful when performing periodic system calculations with variable cells which need to be compared each other. In this way the number of degrees of freedom is kept constant throughout the various simuilations.
- spin_polarize(inp, mpol=1)[source]
Add a collinear spin polarization to the system.
- Parameters:
mpol (int) – spin polarization in Bohr magneton units.
- apply_electric_field(inp, elecfield=[0, 0, 0.001])[source]
Apply an external electric field on the system
- charge_and_polarize(inp)[source]
- Charge the system by removing one electron. Assume that the original
system is closed shell, thus polarize.
- set_scf_method(inp, method='dirmin', mixing_on='density', mixing_scheme='Pulay')[source]
Set the algorithm for scf.
- Parameters:
method (str) – The algoritm chosen. Might be different for the cubic (CS) or linear scaling (LS) algorithm. * dirmin: Direct minimization approach (only CS) * mixing: Mixing scheme (only CS) * hybrid: outer loop of the LS mode * two_levels: two level of accuracy in the scf
mixing_on (str) – May be
"density"
or"potential"
in the"mixing"
case, decide to which quantity the mixing to be performedmixing_scheme (str) –
May be:
Pulay : DIIS mixing on the last 7 iterations
Simple: Simple mixing
Anderson: Anderson scheme
Anderson2: Anderson scheme based on the two pervious iterations
- CG: Conjugate Gradient based on the minimum of the energy with
respect of the potential
Warning
Only the FOE method exhibit asymptotic linear scaling regime.
Todo
Check if the linear scaling case needs another input variable for the mixing of the potential (density)
- add_empty_scf_orbitals(inp, norbs=10)[source]
Insert
norbs
empty orbitals in the scf procedure- Parameters:
norbs (int) – Number of empty orbitals
Warning
In linear scaling case, this is only meaningful for the direct minimization approach.
- write_cubefiles_around_fermi_level(inp, nplot=1)[source]
Writes the
nplot
orbitals around the fermi level in cube format.- Parameters:
nplot (int) – the number of orbitals to print around the fermi level.
Warning
This would work only for the cubic scaling code at present.
- write_orbitals_on_disk(inp, format='binary')[source]
Set the code to write the orbitals on disk in the provided format
- Parameters:
format (str) – The format to write the orbitals with. Accepts the strings: * ‘binary’ * ‘text’ * ‘etsf’ (requires etsf-io enabled)
Todo
Verify if this option works for a linear scaling calulation.
- write_support_functions_on_disk(inp, format='text')[source]
Dump the support functions.
Write the support function which are expressed in wavelet basis in the code as files at the end of the calculation.
- Parameters:
format (str) – the format of the data, can be ‘text’, ‘ETSF’ or ‘binary’ or one of the allowed codes of the output_wf key.
- write_support_function_matrices(inp, format='text')[source]
Write the matrices of the linear scaling formats.
- Parameters:
format (str) –
The format to write the orbitals with. Accepts the strings:
’binary’
’text’
’text_serial’
Todo
Verify if the binary format is available and set the appropriate values
- set_atomic_positions(inp, posinp=None)[source]
Insert the atomic positions as a part of the input dictionary
- read_orbitals_from_disk(inp, mode='cubic')[source]
Read the orbitals from data directory, if available.
- Parameters:
mode (str) – can be ‘cubic’ or ‘linear’. In the first case, orbitals are read as KS objects, whereas in the second kernel (or coeffs) and support functions are expressed.
- set_kernel_guess(inp, mode='kernel')[source]
Method for guessing the kernel at restart.
- Parameters:
mode (str) – Guessing method. Can be: * kernel * coefficients * random * diag * weight * ao * full_kernel * full_coefficients
- set_electronic_temperature(inp, kT=0.001, T=0)[source]
Define the electronic temperature, in AU (
kT
) or K (T
)
- optimize_geometry(inp, method='FIRE', nsteps=50, betax=4.0, frac_fluct=1.0, forcemax=0)[source]
Optimize the geometry of the system
- Parameters:
nsteps (int) – maximum number of atomic steps.
method (str) –
Geometry optimizer. Available keys:
SDCG: A combination of Steepest Descent and Conjugate Gradient
VSSD: Variable Stepsize Steepest Descent method
LBFGS: Limited-memory BFGS
BFGS: Broyden-Fletcher-Goldfarb-Shanno
PBFGS: Same as BFGS with an initial Hessian from a force field
DIIS: Direct inversion of iterative subspace
FIRE: Fast Inertial Relaxation Engine, described by Bitzek et al.
SQNM: Stabilized quasi-Newton minimzer
betax (float) – the step size for the optimization method. This stepsize is system dependent and it has therefore to be determined for each system.
frac_fluct (float) – Fraction of force fluctuations. Stop if fmax < forces_fluct*frac_fluct.
forcemax (float) – Max forces criterion when stop.
- set_xc(inp, xc='PBE')[source]
Set the exchange and correlation approximation
- Parameters:
xc (str) – the Acronym of the XC approximation
Todo
Insert the XC codes corresponding to
libXC
conventions- write_density_on_disk(inp)[source]
Write the charge density on the disk after the last scf convergence
- use_gpu_acceleration(inp, flavour='CUDA')[source]
Employ gpu acceleration when available, for convolutions and Fock operator
- Parameters:
flavour (str) – can be CUDA or OCL. In one case it activates aceleration for the convolutions (useful for dense electronic systems with many k-points). CUDA is an acceleration useful for exact exchange calculations.
- set_scf_iterations(inp, nit=[50, 1])[source]
Set the number of the iteration per loop
- Parameters:
nit (int,list) – integer of the number of iterations. Might be a scalar or a list, up to length two. The first element of the list contains the number of iterations of the direct minimization loop. if
nit
is a scalar, only this contribution is taken into account. The second element is the number of subspace iterations where the hamiltonian is diagonalized in the subspace. For a LS calculation, the number of iteration correspond to the levels of the outer loop.
- change_data_directory(inp, name='')[source]
Modify the name of the
data-
directory. Useful to grab the orbitals from another directory than the run name- Parameters:
name (str) – the name of the run
- calculate_tddft_coupling_matrix(inp, tda=False, rpa=True, fxc=True)[source]
Perform a Casida TDDFT coupling matrix extraction.
- Parameters:
tda (bool) – when
True
, Tamm-Dancoff approximation is used for the extraction of the coupling matrixrpa (bool) – when
False
, the calculation of the RPA term (the linear response of the hartree potential) is switched offfxc (bool) – when
False
, the calculation of the fxc term (the linear response of the XC operator) is switched off.
Note
The arguments
fxc
andrpa
should not be simultaneouslyFalse
.Warning
Presently the LR-TDDFT casida fxc is only available for LDA functionals in ABINIT flavour.
- extract_virtual_states(inp, nvirt=8, davidson=False, norbv=None, itermax_virt=150)[source]
Extract a given number of empty states after the scf cycle.
- connect_run_data(inp, log=None)[source]
Associate the data of the run of a given logfile to the input by retrieving the data directory name of the logfile.
- Parameters:
log (Logfile) – instance of a Logfile class
- calculate_dipole(inp)[source]
Extract the dipole moment from the total charge density.
Note
This function is useful for the linear scaling setup as the cubic scaling approach always calculates the charge density multipoles.
- set_external_potential(inp, mm_pot)[source]
Set the external potential to which the system is submitted outside the QM region
- Parameters:
mm_pot (dict) – dictionary of the external potential which contains the information on the counter-ions
- set_implicit_solvent(inp, itermax=20, minres=0.0001, solvent='water')[source]
- Add an implicit solvent around the system with the
soft-sphere cavity method
- set_dispersion_correction(inp)[source]
Add Grimme’s D3 correction to the energy and the forces
Warning
Only works with Free Boundary conditions
- set_psp_file(inp, filename=None, element=None)[source]
Employ the given PSP file for the provided element
- set_psp_directory(inp, directory='.', elements=None)[source]
Employs all the “psppar.*” files of the given directory as pseudopotentials
- set_psp_nlcc(inp, elements=None)[source]
Employed the built in Non Linear Core Correction (NLCC) pseudopotentials.
- set_psp_krack(inp, functional='PBE', elements=None)[source]
Employed the built in Krack pseudopotentials.
Warning: For certain elements, there are multiple Krack pseudopotentials available with different numbers of electrons. We have choosen the smallest number of electrons as the default behavior. Importantly, the number of electrons when using the LDA or PBE versions may be different.
- optimize_kernel(inp, method='DIAG', dynamic_convergence=True, nit=5, rpnrm=1e-10, alphamix=0.5)[source]
Methods for scf of the density kernel.
Strategies for the optimization of the kernel.
- Parameters:
dynamic_convergence (bool) – if False, the threshold for the convergence of the kernel are not dynamically adjusted
nit (int) – number of scf iterations
rpnrm (float) – convergence criterion, change in density/potential
method (str) – optimization method (‘DIRMIN’, ‘DIAG’, ‘NTPOLY’, ‘FOE’)
alphamix (float) – mixing parameter
- optimize_support_functions(inp, nit=1, gnrm=0.01, diis_history=0)[source]
Methods for scf of the Support functions.
- set_orbital_occupancy(inp, avg={}, up={}, down={}, kpoint=1)[source]
Control the occupation number of the KS orbitals.
With this funtionality one can fix the value of the occupation numbers of the KS orbitals.
- Parameters:
avg (dict) – dictionary of the form {i: f} where i is an integer stating the non-default value of the occupation number, for spin averaged occupancy
up (dict) – same as avg but for the up spin channel
down (dict) – same as avg but for the down spin channel
kpoint (int) – label of the kpoint associated to the occupancy
- set_atomic_occupancy(inp, element=None, iat=None, occupancy={})[source]
Control the occupation number of atomic input guess.
With this funtionality one can fix the value of the occupation numbers of the shells associated to the AO functions employed for the input guess.
- calculate_pdos(inp)[source]
Calculate the Partial Density of states.
Calculate the data needed for Projected DoS on the set of LCAO input wavefunction in the case of the Cubic-Scaling algorithm and onto the Support Functions in the case of the linear scaling.
- calculate_multipoles(inp, yes=True)[source]
Calculate the multipoles on the atoms based on SF.
- Parameters:
yes (bool) – switch off the calculation if False.
- add_cdft_constraint(inp, constraint=None, homo_per_fragment={})[source]
Include a cdft constraint into the SCF cycle.
Useful for Linear scaling formalism.
- Parameters:
constraint (BigDFT.CDFT.CDFTConstraint) – the instantiated constraint
homo_per_fragment (dict) – integer representing the orbital id of the homo of any of the fragments included in the constraint.